Solver API
A Model contains Domains, each of which contain Variables defining Fields, and Reactions with ReactionMethods that operate on the Fields to calculate model time evolution.
Create and initialise
PALEOboxes.create_model_from_config — Methodcreate_model_from_config(
config_file::AbstractString, configmodel::AbstractString;
modelpars::Dict = Dict()
) -> model::Model
create_model_from_config(
config_files,
configmodel::AbstractString;
modelpars::Dict = Dict()
) -> model::ModelConstruct model from a single YAML config_file, or from a collection of config_files, which are read in order and concatenated before being parsed as yaml.
Optional argument modelpars provides parameters that override those in <configmodel>: parameters: section.
PALEOboxes.check_variable_links — Functioncheck_variable_links(model, modeldata; [throw_on_error=true] [, expect_hostdep_varnames=["global.tforce"]]) -> links_ok::BoolCheck all Variables linked correctly, by checking that there are no unexpected host-dependent non-state Variables (ie unlinked Variables)
PALEOboxes.create_modeldata — Functioncreate_modeldata(model::Model [; arrays_eltype=Float64]) -> modeldata::ModelDataCreate a new ModelData struct for model variables with data arrays element type arrays_eltype.
PALEOboxes.ModelData — TypeModelData(model::Model; arrays_eltype::DataType=Float64, allocatenans::Bool=true)Create a ModelData struct containing model data arrays.
One set of data arrays is created with eltype=arrays_eltype, accessed with arrays_idx=1
Additional sets of data arrays may be added by push_arrays_data!, eg in order to support automatic differentiation which requires Dual numbers as the array element type.
Fields
cellranges_all::Vector{AbstractCellRange}: default cellranges covering all domainsdispatchlists_all: default dispatchlists covering all domainssolver_view_all: optional untyped context field for use by external solvers.
PALEOboxes.add_arrays_data! — Functionadd_arrays_data!(
model, modeldata, arrays_eltype::DataType, arrays_tagname::AbstractString;
[method_barrier=nothing] [, generated_dispatch=true] [, kwargs...])Add a data array set to modeldata, allocate memory, and initialize reactiondata.
Element type and tag name are set by arrays_eltype, arrays_tagname
See allocate_variables!(model::Model, modeldata::AbstractModelData, arrays_idx::Int) and [initialize_reactiondata!] for keyword arguments.
PALEOboxes.push_arrays_data! — Functionpush_arrays_data!(modeldata, arrays_eltype::DataType, arrays_tagname::AbstractString)Add an (unallocated) additional array set with element type arrays_eltype.
PALEOboxes.allocate_variables! — Functionallocate_variables!(
vars, modeldata, arrays_idx;
[, eltypemap::Dict{String, DataType}],
[, default_host_dependent_field_data=nothing],
[, allow_base_link=true],
[. use_base_transfer_jacobian=true],
[, use_base_vars=String[]],
[, check_units_opt=:no])Allocate or link memory for VariableDomains vars in modeldata array set arrays_idx
Element type of allocated Arrays is determined by eltype(modeldata, arrays_idx) (the usual case, allowing use of AD types), which can be overridden by Variable :datatype attribute if present (allowing Variables to ignore AD types). :datatype may be either a Julia DataType (eg Float64), or a string to be looked up in eltypemap.
If allow_base_link==true, and any of the following are true a link is made to the base array (arrays_idx=1), instead of allocating a new array in array set arrays_idx: - Variable element type matches modeldata base eltype (arraysidx=1) - `usebasetransferjacobian=trueand Variable:transferjacobianattribute is set - Variable full name is inusebase_vars`
Field data type is determined by Variable :field_data attribute, optionally this can take a default_host_dependent_field_data default for Variables with host_dependent(v)==true (these are Variables with no Target or no Property linked, intended to be external dependencies supplied by the solver).
If check_units_opt != :no then the :units field of linked variable is checked, resulting in either a warning (if check_units_opt=:warn) or error (if check_units_opt=:error).
allocate_variables!(model, modeldata, arrays_idx; kwargs...)Allocate memory for Domain variables for every Domain in model.
See allocate_variables!(domain::Domain, modeldata::AbstractModelData, arrays_idx::Int).
allocate_variables!(domain, modeldata, arrays_idx; [hostdep=false] [, kwargs...])Allocate memory for Domain Variables.
If hostdep=false, only internal Variables are allocated, allowing host-dependent Variables (usually state Variables and derivatives + any external dependencies) to be set to views on host-managed arrays.
See allocate_variables!(vars, modeldata::AbstractModelData, arrays_idx::Int).
PALEOboxes.check_ready — Methodcheck_ready(model, modeldata; [throw_on_error=true]) -> ready::BoolCheck all variable pointers set (ie all arrays allocated for variable data)
PALEOboxes.check_configuration — Methodcheck_configuration(model; [throw_on_error=true]) -> config_ok::BoolCalls optional Reaction check_configuration methods to perform additional configuration checks.
PALEOboxes.initialize_reactiondata! — Methodinitialize_reactiondata!(model::Model, modeldata::AbstractModelData; kwargs...)Processes VarList_...s from ReactionMethods and populates modeldata.sorted_methodsdata_... with sorted lists of ReactionMethods and corresponding Variable accessors.
Optionally calls create_dispatch_methodlists(model, modeldata, modeldata.cellranges_all) to set modeldata.dispatchlists_all to default ReactionMethodDispatchLists for entire model.
Keyword arguments
arrays_indices=1:num_arrays(modeldata):modeldataarrays_idxto generate dispatch lists forcreate_dispatchlists_all=false: true to setmodeldata.dispatchlists_allgenerated_dispatch=true: true to use autogenerated code formodeldata.dispatchlists_all(fast dispatch, slow compile)
PALEOboxes.dispatch_setup — Functiondispatch_setup(model, attribute_name, modeldata, cellranges=modeldata.cellranges_all)Call setup methods, eg to initialize data arrays (including state variables).
attribute_name defines the setup operation performed. dispatch_setup should be called in sequence with attribute_name = :
:setup: initialise Reactions and set up any non-state Variables (eg model grid Variables) (applied tomodeldataarrays_idx=1, values then copied to otherarrays_idx):norm_value: set state Variable values from:norm_valueattribute in .yaml file, and initialise any Reaction state that requires this value (arrays_idx1 only):initial_value(optional): set state Variable values from:initial_valueattribute in .yaml file (arrays_idx1 only)
PALEOboxes.create_dispatch_methodlists — Functioncreate_dispatch_methodlists(model::Model, modeldata::AbstractModelData, arrays_idx::Int, cellranges; kwargs)
-> DispatchMethodLists(list_initialize, list_do)Compile lists of initialize and do methods + corresponding cellrange for main loop do_deriv.
Subset of methods and cellrange to operate on are generated from supplied cellranges.
Keyword arguments
verbose=false: true for additional log outputgenerated_dispatch=true:trueto createReactionMethodDispatchLists (fast dispatch using generated code, slow compile time),falseto createReactionMethodDispatchListNoGen(slow dynamic dispatch, fast compile time)
PALEOboxes.ReactionMethodDispatchList — TypeReactionMethodDispatchListDefines a list of ReactionMethod with corresponding CellRange and views on Variable data (sub)arrays.
Attaching numerical solvers
High-level access to aggregated collections of Variables is provided by VariableAggregator and VariableAggregatorNamed (see Accessing model objects for low-level access).
PALEOboxes.VariableAggregator — TypeVariableAggregator(vars, cellranges, modeldata, arrays_idx) -> VariableAggregatorAggregate multiple VariableDomains into a flattened list (a contiguous Vector).
Creates a VariableAggregator for collection of Variables vars, with indices from corresponding cellranges, for data arrays in modeldata with arrays_idx.
cellranges may contain nothing entries to indicate whole Domain.
This is mostly useful for aggregating state Variables, derivatives, etc to implement an interface to a generic ODE/DAE etc solver.
Values may be copied to and from a Vector using copyto!
PALEOboxes.get_indices — Functionget_indices(va::VariableAggregator, varnamefull::AbstractString; allow_not_found=false) -> indices::UnitRange{Int64}Return indices in flattened Vector corresponding to Variable varnamefull
Base.copyto! — Methodcopyto!(dest::VariableAggregator, src::AbstractVector; sof=1) -> num_copied::IntSet aggregated Variables dest = (contiguous) Vector src.
Optional sof sets first index in src to use.
Base.copyto! — Methodcopyto!(dest::AbstractVector, va::VariableAggregator; dof=1) -> num_copied::IntSet (contiguous) Vector dest = aggregated Variables src
Optional dof sets first index in dest
PALEOboxes.VariableAggregatorNamed — TypeVariableAggregatorNamed(modeldata, arrays_idx=1) -> VariableAggregatorNamed
VariableAggregatorNamed(vars, modeldata, arrays_idx=1) -> VariableAggregatorNamed
VariableAggregatorNamed(va::VariableAggregator; ignore_cellranges=false) -> VariableAggregatorNamedAggregate VariableDomains into nested NamedTuples, with Domain and Variable names as keys and data arrays (from get_data) as values.
Any / characters in Variable names are replaced with __ (double underscore)
This provides direct access to Variables by name, and is mostly useful for testing or for small models.
Fields
vars: nested NamedTuples (domainname, varname) of VariableDomainsvalues: nested NamedTuples (domainname, varname) of data arrays.
Aggregated collections of a subset of Parameters as a flattened Vector (eg for sensitivity studies) is provided by ParameterAggregator:
PALEOboxes.ParameterAggregator — TypeParameterAggregator(parfullnames::Vector{String}, model; eltype=Float64) -> ParameterAggregatorRepresent a subset of model parameters given by parfullnames as a flattened Vector
parfulnames is a Vector of form ["domainname.reactionname.parname", ...] defining a subset of model parameters (NB: must be of type ParDouble or ParDoubleVec ie scalar or vector of Float64).
norm_values can be used to specify normalisation of the flattened parameter vector (defaults to 1.0).
The parameters can then be set from and copied to a flattened Vector using:
copyto!(pa::ParameterAggregator, newvalues::Vector) # set from newvalues .* norm_values
copyto!(currentvalues::Vector, pa::ParameterAggregator) # copy to currentvalues, dividing by norm_values
get_currentvalues(pa::ParameterAggregator) -> currentvalues::VectorThe subset of parameters are then defined by the p parameter Vector used by SciML solvers, and combined with the full set (from the yaml file) to eg solve an ODE to enable sensitivity studies.
eltype can be eg a Dual number to support ForwardDiff automatic differentiation for parameter Jacobians.
Defining CellRanges
PALEOboxes.AbstractCellRange — TypeAbstractCellRangeDefines a range of cells within a Domain.
Fields
All implementations should define:
domain::Domain: theDomaincovered by this cellrange.operatorID::Int: IfoperatorID==0, call allReactions, otherwise only call those with matchingoperatorID(this enables operator splitting).indices: an iterable list of cell indices.
And then may provide subtype-specific fields defining additional ranges of cells.
PALEOboxes.CellRange — TypeCellRange <: AbstractCellRangeDefines a range of cells in a specified Domain as a linear list.
Fields
domainoperatorIDindices: may be any valid Julia indexing range thing eg 1:100, [1 2 3 4], etc
PALEOboxes.CellRangeColumns — TypeCellRangeColumns <: AbstractCellRangeDefines a range of cells in a specified Domain, organised by columns.
Fields
domainoperatorIDindices: iterator through all cells in arbitrary ordercolumns: iterator through columns: columns[n] returns a Pair icol=>cells where cells are ordered top to bottom
PALEOboxes.create_default_cellrange — Functioncreate_default_cellrange(model::Model [; operatorID=0]) -> Vector{AbstractCellRange}Create a Vector of CellRange instances covering the entire model.
PALEOboxes.Grids.create_default_cellrange — Functioncreate_default_cellrange(domain, grid, [; operatorID=0]) -> CellRangeCreate a CellRange for entire domain and supplied operatorID
create_default_cellrange(domain, grid::Nothing [; operatorID=0]) -> CellRangeCreate a CellRange for entire domain. Fallback for a domain with no grid
create_default_cellrange(domain, grid::UnstructuredVectorGrid [; operatorID=0]) -> CellRangeCreate a CellRange for entire domain - use linear index.
create_default_cellrange(domain, grid::UnstructuredColumnGrid [; operatorID=0]) -> CellRangeColumnsCreate a CellRange for entire domain. Return a PB.CellRangeColumns with iterators for columns and cells.
create_default_cellrange(domain, grid::Union{CartesianLinearGrid, CartesianArrayGrid} [; operatorID=0]) -> CellRangeColumnsCreate a CellRange for entire domain. Return a PB.CellRangeColumns provided by cellrange_cartesiantile
PALEOboxes.Grids.cellrange_cartesiantile — Functioncellrange_cartesiantileCreate a range of cells within a region of CartesianLinearGrid specified by rangestuple in a specified PB.Domain
rangestuple is a tuple of Cartesian index ranges eg (1:9, :, :) for a 3D grid.
3D case return a CellRangeColumns
return partitioning of a 3D Domain into n tiles
Main loop
PALEOboxes.do_deriv — Functiondo_deriv(dispatchlists, deltat::Float64=0.0)
do_deriv(dispatchlists, pa::ParameterAggregator, deltat::Float64=0.0)Wrapper function to calculate entire derivative (initialize and do methods) in one call. dispatchlists is from create_dispatch_methodlists.
PALEOboxes.dispatch_methodlist — Functiondispatch_methodlist(dl::ReactionMethodDispatchList, deltat::Float64=0.0)
dispatch_methodlist(dl::ReactionMethodDispatchList, pa::ParameterAggregator, deltat::Float64=0.0)
dispatch_methodlist(dl::ReactionMethodDispatchListNoGen, deltat::Float64=0.0)
dispatch_methodlist(dl::ReactionMethodDispatchListNoGen, pa::ParameterAggregator, deltat::Float64=0.0)Call a list of ReactionMethods.
Implementation
As an optimisation, with dl::ReactionMethodDispatchList uses @generated for Type stability and to avoid dynamic dispatch, instead of iterating over lists.
ReactionMethodDispatchList fields are Tuples hence are fully Typed, the @generated function emits unrolled code with a function call for each Tuple element.
Diagnostics
PALEOboxes.show_methods_setup — Functionshow_methods_setup(model::Model)Display ordered list of Reaction setup methods (registered by add_method_setup!, called by dispatch_setup)
PALEOboxes.show_methods_initialize — Functionshow_methods_initialize(model::Model)Display ordered list of Reaction initialize methods (registered by add_method_initialize!, called by do_deriv at start of each model timestep).
PALEOboxes.show_methods_do — Functionshow_methods_do(model::Model)Display ordered list of Reaction do methods (registered by add_method_do!, called by do_deriv for each model timestep).
PALEOboxes.show_variables — Functionshow_variables(obj, ...) -> TableShow all Variables attached to PALEO object obj
PALEOboxes.show_links — Functionshow_links(vardom::VariableDomain)
show_links(model::Model, varnamefull::AbstractString)
show_links(io::IO, vardom::VariableDomain)
show_links(io::IO, model::Model, varnamefull::AbstractString)Display all VariableReactions linked to this VariableDomain
varnamefull should be of form "<domain name>.<variable name>"
Linked variables are shown as "<domain name>.<reaction name>.<method name>.<local name>"
PALEOboxes.show_parameters — Functionshow_parameters(model) -> DataFrameGet parameters for all reactions in model.
Examples:
Show all model parameters using VS Code table viewer:
julia> vscodedisplay(PB.show_parameters(run.model))Write out all model parameters as csv for import to a spreadsheet:
julia> CSV.write("pars.csv", PB.show_parameters(run.model))show_parameters(react::AbstractReaction) -> DataFrame
show_parameters(classname::AbstractString) -> DataFramelist all parameters for a Reaction react instance or classname
Accessing model objects
Model
PALEOboxes.get_domain — Functionget_domain(model::Model, name::AbstractString; allow_not_found=true) -> Domain or nothing
get_domain(model::Model, domainid) -> DomainGet Domain by name (may be nothing if name not matched) or domainid (range 1:num_domains).
PALEOboxes.get_reaction — Methodget_reaction(model::Model, domainname, reactionname; allow_not_found=false) -> Reaction or nothingGet Reaction by domainname and reaction name
PALEOboxes.set_parameter_value! — Functionset_parameter_value!(model::Model, domainname, reactionname, parname, value)Convenience function to set Parameter value.
PALEOboxes.get_parameter_value — Functionget_parameter_value(model::Model, domainname, reactionname, parname) -> valueConvenience function to get Parameter value.
PALEOboxes.set_variable_attribute! — Functionset_variable_attribute!(model::Model, varnamefull, attributename::Symbol, value)Set varnamefull (of form <domain name>.<var name>) attributename to value.
PALEOboxes.get_variable_attribute — Functionget_variable_attribute(model::Model, varnamefull, attributename::Symbol, missing_value=missing) -> attributevalueGet varnamefull (of form <domain name>.<var name>) attributename.
Domains
PALEOboxes.get_variable — Functionget_variable(obj, varname, ...) -> variableGet variable <: VariableBase by name from PALEO object
PALEOboxes.get_variables — Functionget_variables(varlist::AbstractVarList; flatten=true) -> Vector{VariableReaction}Return VariableReaction in varlist.
If flatten=true, a Vector-of-Vectors is flattened.
get_variables(domain; hostdep=nothing, vfunction=VF_Undefined) -> Vector{VariableDomain}Get domain variables, optionally filtering for subsets based on hostdep and :vfunction attribute
get_variables(domain, filter) -> Vector{VariableDomain}Get subset of domain variables where filter(var) == true.
get_variables(method::AbstractReactionMethod; filterfn = v -> true) -> Vector{VariableReaction}Get VariableReactions from method.varlists as a flat Vector, optionally restricting to those that match filterfn
get_variables(reaction, localname) -> Vector{VariableReaction}
get_variables(reaction; [filterfn=v->true]) -> Vector{VariableReaction}Get matching Variables from all ReactionMethods.
PALEOboxes.get_host_variables — Functionget_host_variables(domain, vfunction; [match_deriv_suffix=""] [, operatorID=0] [, exclude_var_nameroots=[]] [, verbose=false])
-> (host_vars, host_deriv_vars)Get state Variables with VariableFunction vfunction, and optionally corresponding time derivative with VariableFunction VF_Deriv and name matching hostvarname*<match_deriv_suffix`>.
Optionally filter by operatorID, omit Variables with name matching exclude_var_nameroots.
PALEOboxes.get_unallocated_variables — Functionget_unallocated_variables(domain, modeldata, arrays_idx::Int) -> Vector{VariableDomain}Return any unallocated variables (host-dependent variables which have no data pointer set)
PALEOboxes.get_reaction — Methodget_reaction(domain, reactname; allow_not_found) -> Reaction or nothingGet a reaction by name.
VariableDomain
Low-level access to individual Variables.
PALEOboxes.set_data! — Functionset_data!(var::VariableDomain, modeldata, arrays_idx::Int, data)Set VariableDomain to a newly-created Field containing data.
PALEOboxes.get_data — Functionget_data(var::VariableDomain, modeldata::AbstractModelData, arrays_idx::Int=1) -> field.values
get_data(var::VariableDomain, domaindata::AbstractDomainData) -> field.valuesGet Variable var data array from Field.values
get_data(va::VariableAggregator) -> VectorAllocate Vector and set to values of aggregated Variables va.